linalg_svd Module



Interfaces

public interface svd

  • private module subroutine svd_dbl(a, s, u, vt, work, olwork, err)

    Computes the singular value decomposition of an M-by-N matrix such that where is an M-by-M orthogonal matrix, is an M-by-N diagonal matrix containing the singular values, and is an N-by-N orthogonal matrix.

    Arguments

    Type IntentOptional Attributes Name
    real(kind=real64), intent(inout), dimension(:,:) :: a

    On input, the M-by-N matrix to factor. The matrix is overwritten on output.

    real(kind=real64), intent(out), dimension(:) :: s

    A MIN(M, N)-element array containing the singular values of a sorted in descending order.

    real(kind=real64), intent(out), optional, dimension(:,:) :: u

    An optional argument, that if supplied, is used to contain the orthogonal matrix from the decomposition. The matrix contains the left singular vectors, and can be either M-by-M (all left singular vectors are computed), or M-by-MIN(M,N) (only the first MIN(M, N) left singular vectors are computed).

    real(kind=real64), intent(out), optional, dimension(:,:) :: vt

    An optional argument, that if supplied, is used to contain the transpose of the N-by-N orthogonal matrix . The matrix contains the right singular vectors.

    real(kind=real64), intent(out), optional, target, dimension(:) :: work

    An optional input, that if provided, prevents any local memory allocation. If not provided, the memory required is allocated within. If provided, the length of the array must be at least olwork.

    integer(kind=int32), intent(out), optional :: olwork

    An optional output used to determine workspace size. If supplied, the routine determines the optimal size for work, and returns without performing any actual calculations.

    class(errors), intent(inout), optional, target :: err

    The error object to be updated.

  • private subroutine svd_cmplx(a, s, u, vt, work, olwork, rwork, err)

    Computes the singular value decomposition of an M-by-N matrix such that where is an M-by-M orthogonal matrix, is an M-by-N diagonal matrix containing the singular values, and is an N-by-N orthogonal matrix.

    Arguments

    Type IntentOptional Attributes Name
    complex(kind=real64), intent(inout), dimension(:,:) :: a

    On input, the M-by-N matrix to factor. The matrix is overwritten on output.

    real(kind=real64), intent(out), dimension(:) :: s

    A MIN(M, N)-element array containing the singular values of a sorted in descending order.

    complex(kind=real64), intent(out), optional, dimension(:,:) :: u

    An optional argument, that if supplied, is used to contain the orthogonal matrix from the decomposition. The matrix contains the left singular vectors, and can be either M-by-M (all left singular vectors are computed), or M-by-MIN(M,N) (only the first MIN(M, N) left singular vectors are computed).

    complex(kind=real64), intent(out), optional, dimension(:,:) :: vt

    An optional argument, that if supplied, is used to contain the conjugate transpose of the N-by-N orthogonal matrix . The matrix contains the right singular vectors.

    complex(kind=real64), intent(out), optional, target, dimension(:) :: work

    An optional input, that if provided, prevents any local memory allocation. If not provided, the memory required is allocated within. If provided, the length of the array must be at least olwork.

    integer(kind=int32), intent(out), optional :: olwork

    An optional output used to determine workspace size. If supplied, the routine determines the optimal size for work, and returns without performing any actual calculations.

    real(kind=real64), intent(out), optional, target, dimension(:) :: rwork

    An optional input, that if provided, prevents any local memory allocation for real-valued workspaces. If not provided, the memory required is allocated within. If provided, the length of the array must be at least 5 * MIN(M, N).

    class(errors), intent(inout), optional, target :: err

    The error object to be updated.